In [39]:
# necessary packages #

#using Pkg
#Pkg.add("Distances")
using Distributions
using Random
using Distances
using LinearAlgebra
using SparseArrays
using IterativeSolvers
using ProgressMeter
In [40]:
include("util.j")
Out[40]:
getAD (generic function with 1 method)
In [77]:
# unnecessary packages #

#using Pkg
#Pkg.add("UnicodePlots")
using UnicodePlots   # check the structure of the sparse matrix
using BenchmarkTools

using StatsPlots
using MCMCChains
using PrettyTables
In [42]:
#using Pkg
#Pkg.add("ProgressMeter");
In [43]:
# Set the parameters for SLMC model #

N = 1200 # sample size
N1 = 1000; N2 = 1000;
q = 2; p = 2; K = 3
Σ = [0.3 0.1
     0.1 0.2];
β = [1.0 -1.0
     -5.0 2.0];
ϕ1 = 5.0; ϕ2 = 12.0; ϕ3 = 30.0; ν1 = 0.5; ν2 = 0.5; ν3 = 0.5; # parameter for the independent F
Λ = [3.0 0.0
     -2.0 3.0
     1.0 -2.0]; # loading matrix

# priors #
μβ = fill(0.0, p, q);  =[[100.0 0.0]; [0.0 100.0]];
μΛ = fill(0.0, K, q);  =[[100.0 0.0 0.0]; [0.0 100.0 0.0]; [0.0 0.0 100.0]];
νΣ = 2 * q; ΨΣ = [[1.0 0.0]; [0.0 1.0]];
ϕU = 300 / sqrt(2); ϕL = 3 / sqrt(2);
In [44]:
# Generate simulation data #

Random.seed!(1234);
coords = rand(2, N);                                          # random location over unit square
X = hcat(fill(1, (N,)), rand(N));                             # design matrix
D = pairwise(Euclidean(), coords, dims = 2);                  # distance matrix
ρ1 = exp.(-ϕ1 * D); ρ2 = exp.(-ϕ2 * D);ρ3 = exp.(-ϕ3 * D);    # covariance matrix
ω = [rand(MvNormal(ρ1), 1) rand(MvNormal(ρ2), 1) rand(MvNormal(ρ3), 1)] * Λ; # latent process
Y = X * β + ω + transpose(rand(MvNormal(Σ), N));              # response matrix
In [45]:
# Some data preparations #

ordx = sortperm(coords[1, :]);                                # sort order based on the first coordinates
X_ord = X[ordx, :]; Y_ord = Y[ordx, :]; ω_ord = ω[ordx, :];   # sorted data
ω_incp_obs = ω_ord + fill(1.0, (N, 1)) * transpose(β[1, :]); # latent process + intercept
coords_ord = coords[:, ordx];
S1_ind = sample(1:N, N1, replace = false, ordered = true);    # observed location index for 1st response
S2_ind = sample(1:N, N2, replace = false, ordered = true);    # observed location index for 2nd response
S = sort(union(S1_ind, S2_ind));                              # observed index set
M1_ind = setdiff(S, S1_ind);                                  # in S not in S1
M2_ind = setdiff(S, S2_ind);                                  # in S not in S2 
obs_ind = vcat(S1_ind, S2_ind .+ N);              # index of the observed location for all response among N locations
perm_ind = sortperm(vcat(S1_ind, S2_ind));                    # the vector of the permutation 

v1 = zeros(N); v1[S1_ind] .= 1;
v2 = zeros(N); v2[S2_ind] .= 1;
index_S = (2^0 * v1 + 2^1 * v2);                              # build index indicating which response are observed
M1_Sind = findall(x -> x == 2^1, index_S[S]);                 # index of M1 among S
M2_Sind = findall(x -> x == 2^0, index_S[S]);                 # index of M2 among S

m = 10; n = length(S);                                        # number of nearest neighbor                       
NN = BuildNN(coords_ord[:, S], m, 1.0);                            # build nearest neighbor 
nnIndx_col = vcat(NN.nnIndx, 1:n);                            # the index of columns
nnIndx_row = zeros(Int64, 0);                                               
for i in 2:m
    nnIndx_row = vcat(nnIndx_row, fill(i, i-1))
end
nnIndx_row = vcat(nnIndx_row, repeat((m + 1):n, inner = m), 1:n);  # the index of rows

 = cholesky();  = cholesky(); 
RWM_scale = 0.2;                                              # random-walk metropolis step size scale
In [46]:
# check the plot of the data 
using RCall
@rput ω_incp_obs
@rput coords_ord
@rput S
R"""
library(MBA)
library(classInt)
library(RColorBrewer)
library(sp)
library(coda)
library(spBayes)
library(fields)

h <- 12
surf.raw1 <- mba.surf(cbind(t(coords_ord[, S]), ω_incp_obs[S, 1]), no.X = 300, no.Y = 300, 
              exten = TRUE, sp = TRUE, h = h)$xyz.est
surf.raw2 <- mba.surf(cbind(t(coords_ord[, S]), ω_incp_obs[S, 2]), no.X = 300, no.Y = 300, 
              exten = TRUE, sp = TRUE, h = h)$xyz.est
surf.brks <- classIntervals(c(surf.raw1[["z"]], surf.raw2[["z"]]), 500, 'pretty')$brks
col.pal <- colorRampPalette(brewer.pal(11,'RdBu')[11:1])
xlim <- c(0, 1.13)

zlim <- range(c(surf.raw1[["z"]], surf.raw2[["z"]]))

# size for the mapping of w               
width <- 360
height <- 360
pointsize <- 16

png(paste("./pics/map-w1_incp-true.png", sep = ""), 
    width = width, height = height, pointsize = pointsize, family = "Courier")
par(mfrow = c(1, 1))
##Obs
i <- as.image.SpatialGridDataFrame(surf.raw1)
plot(t(coords_ord[, S]), typ="n", cex=0.5, xlim=xlim, axes=FALSE, ylab="y", xlab="x") 
     #main = "true")
axis(2, las=1)
axis(1)
image.plot(i, add=TRUE, col=rev(col.pal(length(surf.brks)-1)), zlim=zlim)
dev.off()

png(paste("./pics/map-w2_incp-true.png", sep = ""), 
    width = width, height = height, pointsize = pointsize, family = "Courier")
par(mfrow = c(1, 1))
##Obs
i <- as.image.SpatialGridDataFrame(surf.raw2)
plot(t(coords_ord[, S]), typ="n", cex=0.5, xlim=xlim, axes=FALSE, ylab="y", xlab="x") 
     #main = "true")
axis(2, las=1)
axis(1)
image.plot(i, add=TRUE, col=rev(col.pal(length(surf.brks)-1)), zlim=zlim)
dev.off()

"""
Out[46]:
RObject{IntSxp}
null device 
          1 
In [47]:
# preallocation #

#F = Array{Float64,2}(undef, n , 3);                           # preallocate the matrix F
μ_m1 = Array{Float64, 2}(undef, length(M1_ind), q);
μ_m2 = Array{Float64, 2}(undef, length(M2_ind), q);
nIndx = length(NN.nnIndx);
A1 = Array{Float64}(undef, nIndx); D1 = Array{Float64}(undef, n);
A2 = Array{Float64}(undef, nIndx); D2 = Array{Float64}(undef, n);
A3 = Array{Float64}(undef, nIndx); D3 = Array{Float64}(undef, n);
I_A1 = SparseMatrixCSC{Float64,Int64};
I_A2 = SparseMatrixCSC{Float64,Int64}; 
I_A3 = SparseMatrixCSC{Float64,Int64}; 
A1_new = Array{Float64}(undef, nIndx); D1_new = Array{Float64}(undef, n);
A2_new = Array{Float64}(undef, nIndx); D2_new = Array{Float64}(undef, n);
A3_new = Array{Float64}(undef, nIndx); D3_new = Array{Float64}(undef, n);
I_A1_new = SparseMatrixCSC{Float64,Int64};
I_A2_new = SparseMatrixCSC{Float64,Int64}; 
I_A3_new = SparseMatrixCSC{Float64,Int64}; 
Ystar = vcat(Y_ord[S, :], .L \ μβ, .L \ μΛ);             # will be updated after imputing missing response
Xstar = vcat([X_ord[S, :] fill(0.0, n, K)], [inv(.L) fill(0.0, p, K)], 
    [fill(0.0, K, p) inv(.L)]);      
Ψstar = fill(0.0, q, q); νstar = νΣ + n;
μγstar = vcat(μβ, μΛ); Vγstar = fill(0.0, p + K, p + K);
Y_Xm = fill(0.0, n, q); invVγstar = fill(0.0, p + K, p + K);

MCMC sampling algorithm Q1: priors for $\nu_i$ Q2: $\phi_i$ may not be consistant, since the order can change

In [48]:
# Preallocation for MCMC samples and Initalization #
N_sam = 20000;
γ_sam = Array{Float64, 3}(undef, p + K, q, N_sam + 1);
Σ_sam = Array{Float64, 3}(undef, q, q, N_sam + 1);
F_sam = Array{Float64, 3}(undef, n, K, N_sam);
Y_m1_sam = Array{Float64, 2}(undef, length(M1_ind), N_sam);
Y_m2_sam = Array{Float64, 2}(undef, length(M2_ind), N_sam);
A_sam = Array{Float64, 1}(undef, N_sam); 
lh_old = 1; lh_new = 1;     # record the likelihood for updating ranges

ϕ_sam = Array{Float64, 2}(undef, K, N_sam + 1);

γ_sam[:, :, 1] = vcat([[0.0 0.0]; [0.0 0.0]], [[1.0 0.0]; [0.0 1.0]; [1.0 1.0]]);
Σ_sam[:, :, 1] = [[0.5 0.1]; [0.1 0.5]];
ϕ_sam[:, 1] = [6, 20, 10];
In [49]:
# for loop for MCMC chain #
prog = Progress(N_sam, 1, "Computing initial pass...", 50)
for l in 1:N_sam
    # Build the matrix D_Sigma_o^{1/2} #
    Dic_diag = Dict(2^0 => sparse(1I, 1, 1) * (1 / sqrt(Σ_sam[:, :, l][1, 1])), 
        2^1 => sparse(1I, 1, 1) * (1 / sqrt(Σ_sam[:, :, l][2, 2])), 
        (2^0 + 2^1)=> sparse(sqrt(inv(Σ_sam[:, :, l]))));
    invD = blockdiag([Dic_diag[i] for i in index_S if i > 0]...);
                    
    # Build the matrix for the first iteration #
    if l == 1
        getAD(coords_ord[:, S], NN.nnIndx, NN.nnDist, NN.nnIndxLU, ϕ_sam[1, l], 0.5, A1, D1);
        getAD(coords_ord[:, S], NN.nnIndx, NN.nnDist, NN.nnIndxLU, ϕ_sam[2, l], 0.5, A2, D2);
        getAD(coords_ord[:, S], NN.nnIndx, NN.nnDist, NN.nnIndxLU, ϕ_sam[3, l], 0.5, A3, D3);
        I_A1 = sparse(nnIndx_row, nnIndx_col, vcat(-A1, ones(n)));
        I_A2 = sparse(nnIndx_row, nnIndx_col, vcat(-A2, ones(n)));
        I_A3 = sparse(nnIndx_row, nnIndx_col, vcat(-A3, ones(n)));
    end
                    

    # Build Ytilde Xtilde
    Ytilde = vcat(invD * vcat(Y_ord[S1_ind, 1] - X_ord[S1_ind, :] * γ_sam[1:p, 1, l], 
                            Y_ord[S2_ind, 2] - X_ord[S2_ind, :] * γ_sam[1:p, 2, l])[perm_ind], zeros(K * n));
    Xtilde = vcat(invD * kron(sparse(transpose(γ_sam[(p + 1):(p + K), :, l])), 
                            sparse(1:N, 1:N, ones(N)))[obs_ind, 
                            vcat(S, S .+ N, S .+ (2 * N))][perm_ind, :],
             blockdiag(Diagonal(1 ./ sqrt.(D1)) * I_A1, Diagonal(1 ./ sqrt.(D2)) * I_A2,
                            Diagonal(1 ./ sqrt.(D3)) * I_A3));
                
    # use LSMR to generate sample of F #       
    nsam = length(Ytilde);
    F_sam[:, :, l] = reshape(lsmr(Xtilde, Ytilde + rand(Normal(), nsam)), :, K);
                    
    # impute missing response  over S#
    Xstar[1:n, (p + 1):(p + K)] = F_sam[:, :, l];        # update matrix Xstar with F
    mul!(μ_m1, Xstar[M1_Sind, :], γ_sam[:, :, l]);
    mul!(μ_m2, Xstar[M2_Sind, :], γ_sam[:, :, l]);

    Y_m1_sam[:, l] = μ_m1[:, 1] + (Σ_sam[1, 2, l] / Σ_sam[2, 2, l]) * 
            (Y_ord[M1_ind, 2] - μ_m1[:, 2]) + 
            rand(Normal(0, sqrt(Σ_sam[1, 1, l] - Σ_sam[1, 2, l]^2 / Σ_sam[2, 2, l])), length(M1_ind));
    Y_m2_sam[:, l] = μ_m2[:, 2] + (Σ_sam[2, 1, l] / Σ_sam[1, 1, l]) * 
            (Y_ord[M2_ind, 1] - μ_m2[:, 1]) + 
            rand(Normal(0, sqrt(Σ_sam[2, 2, l] - Σ_sam[2, 1, l]^2 / Σ_sam[1, 1, l])), length(M2_ind)); # improve?...
    
                    
    # use MNIW to sample γ Σ #
    Ystar[M1_Sind, 1] = Y_m1_sam[:, l];              # update Ystar with imputed response
    Ystar[M2_Sind, 2] = Y_m2_sam[:, l]; 
    invVγstar = cholesky(Xstar'Xstar);
    mul!(μγstar, transpose(Xstar), Ystar); μγstar = invVγstar.U \ (invVγstar.L \ μγstar);
    Y_Xm = BLAS.gemm('N', 'N', -1.0, Xstar, μγstar) + Ystar;
    mul!(Ψstar, transpose(Y_Xm), Y_Xm); BLAS.axpy!(1.0, ΨΣ, Ψstar);

    Σ_sam[:, :, l + 1] = rand(InverseWishart(νstar, Ψstar), 1)[1];    # sample Σ
    γ_sam[:, :, l + 1] = (invVγstar.U \ reshape(rand(Normal(), (p + K) * q), (p + K), q)) * 
                    cholesky(Σ_sam[:, :, l + 1]).U + μγstar;          # sample γ    
                    
                    
    # use metropolis-hasting to update range
    ϕ_sam[:, l + 1] = ϕ_sam[:, l] + RWM_scale * rand(Normal(), K); # propose next sample point
    if all(x -> (x < ϕU && x > ϕL), ϕ_sam[:, l + 1])
        lh_old = -0.5 * (sum(log.(D1)) + sum(log.(D2)) + sum(log.(D3)) + 
                    norm((I_A1 * F_sam[:, 1, l]) ./ sqrt.(D1))^2 + 
                    norm((I_A2 * F_sam[:, 2, l]) ./ sqrt.(D2))^2 + 
                    norm((I_A3 * F_sam[:, 3, l]) ./ sqrt.(D3))^2);
        getAD(coords_ord[:, S], NN.nnIndx, NN.nnDist, NN.nnIndxLU, ϕ_sam[1, l + 1], 0.5, A1_new, D1_new);
        getAD(coords_ord[:, S], NN.nnIndx, NN.nnDist, NN.nnIndxLU, ϕ_sam[2, l + 1], 0.5, A2_new, D2_new);
        getAD(coords_ord[:, S], NN.nnIndx, NN.nnDist, NN.nnIndxLU, ϕ_sam[3, l + 1], 0.5, A3_new, D3_new);
        I_A1_new = sparse(nnIndx_row, nnIndx_col, vcat(-A1_new, ones(n)));
        I_A2_new = sparse(nnIndx_row, nnIndx_col, vcat(-A2_new, ones(n)));
        I_A3_new = sparse(nnIndx_row, nnIndx_col, vcat(-A3_new, ones(n)));            
        lh_new = -0.5 * (sum(log.(D1_new)) + sum(log.(D2_new)) + sum(log.(D3_new)) + 
                    norm((I_A1_new * F_sam[:, 1, l]) ./ sqrt.(D1_new))^2 + 
                    norm((I_A2_new * F_sam[:, 2, l]) ./ sqrt.(D2_new))^2 + 
                    norm((I_A3_new * F_sam[:, 3, l]) ./ sqrt.(D3_new))^2);
        A_sam[l] = exp(lh_new - lh_old);
        if rand(1)[1] < A_sam[l]
            I_A1 = copy(I_A1_new); I_A2 = copy(I_A2_new); I_A3 = copy(I_A3_new); # update and update the corresponding I_A D
            D1 = copy(D1_new); D2 = copy(D2_new); D3 = copy(D3_new); 
        else
            ϕ_sam[:, l + 1] = ϕ_sam[:, l]; # Don't update
        end         
    else
        A_sam[l] = 0.0;
        ϕ_sam[:, l + 1] = ϕ_sam[:, l];   # Don't update when falling out of the supports
    end                      
    
    next!(prog) # monitor the progress
end
Computing initial pass...100%|██████████████████████████████████████████████████| Time: 1:20:57:57

MCMC Chain check

In [50]:
β_pos_sam = Array{Float64, 3}(undef, N_sam + 1, p * q, 1);
β_pos_sam[:, :, 1] = hcat(γ_sam[1, 1, :], γ_sam[1, 2, :], γ_sam[2, 1, :], γ_sam[2, 2, :]);
β_chain = Chains(β_pos_sam);
 = plot(β_chain)
Out[50]:
0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 0 1 2 3 Param1 Iteration Sample value 0 1 2 3 0.0 0.2 0.4 0.6 Param1 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -2.5 -2.0 -1.5 -1.0 -0.5 0.0 Param2 Iteration Sample value -3 -2 -1 0 0.0 0.2 0.4 0.6 0.8 Param2 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -5 -4 -3 -2 -1 0 Param3 Iteration Sample value -5 -4 -3 -2 -1 0 0.0 0.5 1.0 1.5 2.0 2.5 Param3 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 0.0 0.5 1.0 1.5 2.0 2.5 Param4 Iteration Sample value 0.0 0.5 1.0 1.5 2.0 2.5 0.0 0.5 1.0 1.5 2.0 Param4 Sample value Density
In [51]:
Λ_pos_sam = Array{Float64, 3}(undef, N_sam + 1, K * q, 1);
Λ_pos_sam[:, :, 1] = hcat(γ_sam[3, 1, :], γ_sam[3, 2, :], γ_sam[4, 1, :], γ_sam[4, 2, :], 
    γ_sam[5, 1, :], γ_sam[5, 2, :]);
Λ_chain = Chains(Λ_pos_sam);
 = plot(Λ_chain)
Out[51]:
0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 1 2 3 4 Param1 Iteration Sample value 1 2 3 4 5 0.0 0.1 0.2 0.3 0.4 Param1 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -1.5 -1.0 -0.5 0.0 0.5 1.0 1.5 Param2 Iteration Sample value -2 -1 0 1 2 0.0 0.1 0.2 0.3 0.4 0.5 0.6 Param2 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -2.5 -2.0 -1.5 -1.0 -0.5 0.0 Param3 Iteration Sample value -2.5 -2.0 -1.5 -1.0 -0.5 0.0 0.00 0.25 0.50 0.75 1.00 Param3 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 1.0 1.5 2.0 2.5 3.0 3.5 Param4 Iteration Sample value 1.0 1.5 2.0 2.5 3.0 3.5 0.00 0.25 0.50 0.75 1.00 Param4 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -4 -3 -2 -1 0 1 Param5 Iteration Sample value -5 -4 -3 -2 -1 0 1 0.0 0.1 0.2 0.3 0.4 Param5 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 1 2 3 4 Param6 Iteration Sample value 1 2 3 4 0.0 0.1 0.2 0.3 0.4 0.5 0.6 Param6 Sample value Density
In [52]:
ϕ_pos_sam = Array{Float64, 3}(undef, N_sam + 1, K, 1);
ϕ_pos_sam[:, :, 1] = hcat(ϕ_sam[1, :], ϕ_sam[2, :], ϕ_sam[3, :]);
ϕ_chain = Chains(ϕ_pos_sam);
 = plot(ϕ_chain)
Out[52]:
0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 2.5 5.0 7.5 10.0 12.5 Param1 Iteration Sample value 2.5 5.0 7.5 10.0 12.5 0.00 0.05 0.10 0.15 0.20 Param1 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 15 20 25 30 35 40 45 Param2 Iteration Sample value 10 20 30 40 50 0.00 0.02 0.04 0.06 Param2 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 5 10 15 20 25 30 Param3 Iteration Sample value 0 10 20 30 0.000 0.025 0.050 0.075 0.100 Param3 Sample value Density
In [53]:
Σ_pos_sam = Array{Float64, 3}(undef, N_sam + 1, q * q, 1);
Σ_pos_sam[:, :, 1] = hcat(Σ_sam[1, 1, :], Σ_sam[1, 2, :], Σ_sam[2, 1, :], Σ_sam[2, 2, :]);
Σ_chain = Chains(Σ_pos_sam);
 = plot(Σ_chain)
Out[53]:
0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 0.5 1.0 1.5 2.0 Param1 Iteration Sample value 0.0 0.5 1.0 1.5 2.0 0 1 2 3 4 5 Param1 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -1.5 -1.0 -0.5 0.0 Param2 Iteration Sample value -1.5 -1.0 -0.5 0.0 0 2 4 6 8 Param2 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -1.5 -1.0 -0.5 0.0 Param3 Iteration Sample value -1.5 -1.0 -0.5 0.0 0 2 4 6 8 Param3 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 0.0 0.5 1.0 1.5 2.0 Param4 Iteration Sample value 0.0 0.5 1.0 1.5 2.0 0 1 2 3 4 5 6 Param4 Sample value Density
In [78]:
ω_incp_obs_pos_sam = Array{Float64, 3}(undef, n, q, N_sam);
lll = fill(1.0, (n, 1));
for i in 1:N_sam
    ω_incp_obs_pos_sam[:, :, i] = F_sam[:, :, i] * γ_sam[(p + 1):(p + K), :, i + 1] + 
                     lll * transpose(γ_sam[1, :, i + 1]);
end
truncindex = 1;#Integer(trunc(N_sam / 2));
ω_incp_pos_sam = Array{Float64, 3}(undef, N_sam - truncindex  + 1, 3, 1);
ω_incp_pos_sam[:, :, 1] = hcat(ω_incp_obs_pos_sam[1, 1, truncindex:N_sam], 
    ω_incp_obs_pos_sam[1, 2, truncindex:N_sam], ω_incp_obs_pos_sam[200, 1, truncindex:N_sam]);
ω_incp_chain = Chains(ω_incp_pos_sam);
 = plot(ω_incp_chain)
Out[78]:
0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -5 -4 -3 -2 -1 Param1 Iteration Sample value -5 -4 -3 -2 -1 0.0 0.2 0.4 0.6 0.8 Param1 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -3 -2 -1 0 1 Param2 Iteration Sample value -3 -2 -1 0 1 2 0.00 0.25 0.50 0.75 1.00 Param2 Sample value Density 0 5.0×10 3 1.0×10 4 1.5×10 4 2.0×10 4 -1 0 1 2 3 Param3 Iteration Sample value -1 0 1 2 3 0.0 0.2 0.4 0.6 0.8 Param3 Sample value Density

Posterior Inference

In [79]:
N_Inf_burn = Integer(trunc(N_sam / 2));
ω_incp_obs_pos_qt = Array{Float64, 3}(undef, n, q, 3);
for j in 1:q
    for i in 1:n
        ω_incp_obs_pos_qt[i, j, :] = quantile(ω_incp_obs_pos_sam[i, j, N_Inf_burn:N_sam], [0.025, 0.5, 0.975])
    end
end
# count the covarage of 95% CI #
count = fill(0.0, 2);
for j in 1:q
    for i in 1:n
        count[j] = count[j] + 
        ((ω_incp_obs_pos_qt[i, j, 1] < ω_incp_obs[S[i], j]) && 
            (ω_incp_obs_pos_qt[i, j, 3] > ω_incp_obs[S[i], j]))
    end
end
count
Out[79]:
2-element Array{Float64,1}:
 1113.0
 1083.0
In [80]:
count ./ n
Out[80]:
2-element Array{Float64,1}:
 0.9561855670103093
 0.9304123711340206
In [82]:
summary_table = Array{Float64, 2}(undef, (p - 1) * q + q * q, 5);
summary_table[1, :] = vcat(β[2, 1], mean(γ_sam[2, 1, N_Inf_burn:(N_sam + 1)]),
    quantile(γ_sam[2, 1, N_Inf_burn:(N_sam + 1)], [0.5, 0.025, 0.975]));
summary_table[2, :] = vcat(β[2, 2], mean(γ_sam[2, 2, N_Inf_burn:(N_sam + 1)]),
    quantile(γ_sam[2, 2, N_Inf_burn:(N_sam + 1)], [0.5, 0.025, 0.975]));
summary_table[3, :] = vcat(Σ[1, 1], mean(Σ_sam[1, 1, N_Inf_burn:(N_sam + 1)]),
    quantile(Σ_sam[1, 1, N_Inf_burn:(N_sam + 1)], [0.5, 0.025, 0.975]));
summary_table[4, :] = vcat(Σ[1, 2], mean(Σ_sam[1, 2, N_Inf_burn:(N_sam + 1)]),
    quantile(Σ_sam[1, 2, N_Inf_burn:(N_sam + 1)], [0.5, 0.025, 0.975]));
summary_table[5, :] = vcat(Σ[2, 1], mean(Σ_sam[2, 1, N_Inf_burn:(N_sam + 1)]),
    quantile(Σ_sam[2, 1, N_Inf_burn:(N_sam + 1)], [0.5, 0.025, 0.975]));
summary_table[6, :] = vcat(Σ[2, 2], mean(Σ_sam[2, 2, N_Inf_burn:(N_sam + 1)]),
    quantile(Σ_sam[2, 2, N_Inf_burn:(N_sam + 1)], [0.5, 0.025, 0.975]));
summary_table = round.(summary_table; digits = 3 );
rnames = ["β[2, 1]", "β[2, 2]", "Σ[1, 1]", "Σ[1, 2]", "Σ[2, 1]", "Σ[2, 2]"];
summary_table = [rnames summary_table];
pretty_table(summary_table,  ["" "true" "mean" "median" "2.5%" "97.5%"], markdown)
|         | true |   mean | median |   2.5% |  97.5% |
|---------|------|--------|--------|--------|--------|
| β[2, 1] | -5.0 | -4.839 | -4.842 | -5.123 | -4.545 |
| β[2, 2] |  2.0 |  1.683 |  1.682 |  1.343 |  2.018 |
| Σ[1, 1] |  0.3 |  0.283 |  0.279 |  0.151 |  0.443 |
| Σ[1, 2] |  0.1 |  0.015 |  0.019 |  -0.09 |  0.104 |
| Σ[2, 1] |  0.1 |  0.015 |  0.019 |  -0.09 |  0.104 |
| Σ[2, 2] |  0.2 |  0.153 |  0.143 |  0.058 |  0.298 |

Plot the latent process + intercept

In [55]:
# check the plot of the data 
using RCall
@rput ω_incp_obs
@rput coords_ord
@rput S
@rput ω_incp_obs_pos_qt
R"""
library(MBA)
library(classInt)
library(RColorBrewer)
library(sp)
library(coda)
library(spBayes)
library(fields)

h <- 12
surf.raw1 <- mba.surf(cbind(t(coords_ord[, S]), ω_incp_obs[S, 1]), no.X = 300, no.Y = 300, 
              exten = TRUE, sp = TRUE, h = h)$xyz.est
surf.raw2 <- mba.surf(cbind(t(coords_ord[, S]), ω_incp_obs[S, 2]), no.X = 300, no.Y = 300, 
              exten = TRUE, sp = TRUE, h = h)$xyz.est

surf.pos1 <- mba.surf(cbind(t(coords_ord[, S]), ω_incp_obs_pos_qt[, 1, 2]), no.X = 300, no.Y = 300, 
              exten = TRUE, sp = TRUE, h = h)$xyz.est
surf.pos2 <- mba.surf(cbind(t(coords_ord[, S]), ω_incp_obs_pos_qt[, 2, 2]), no.X = 300, no.Y = 300, 
              exten = TRUE, sp = TRUE, h = h)$xyz.est

surf.brks <- classIntervals(c(surf.raw1[["z"]], surf.raw2[["z"]]), 500, 'pretty')$brks
col.pal <- colorRampPalette(brewer.pal(11,'RdBu')[11:1])
xlim <- c(0, 1.13)

zlim <- range(c(surf.raw1[["z"]], surf.raw2[["z"]], surf.pos1[["z"]], surf.pos2[["z"]]))

# size for the mapping of w               
width <- 360
height <- 360
pointsize <- 16

png(paste("./pics/map-w1_incp-true.png", sep = ""), 
    width = width, height = height, pointsize = pointsize, family = "Courier")
par(mfrow = c(1, 1))
##Obs
i <- as.image.SpatialGridDataFrame(surf.raw1)
plot(t(coords_ord[, S]), typ="n", cex=0.5, xlim=xlim, axes=FALSE, ylab="y", xlab="x") 
     #main = "true")
axis(2, las=1)
axis(1)
image.plot(i, add=TRUE, col=rev(col.pal(length(surf.brks)-1)), zlim=zlim)
dev.off()

png(paste("./pics/map-w2_incp-true.png", sep = ""), 
    width = width, height = height, pointsize = pointsize, family = "Courier")
par(mfrow = c(1, 1))
##Obs
i <- as.image.SpatialGridDataFrame(surf.raw2)
plot(t(coords_ord[, S]), typ="n", cex=0.5, xlim=xlim, axes=FALSE, ylab="y", xlab="x") 
     #main = "true")
axis(2, las=1)
axis(1)
image.plot(i, add=TRUE, col=rev(col.pal(length(surf.brks)-1)), zlim=zlim)
dev.off()

png(paste("./pics/map-w1_incp-posm.png", sep = ""), 
    width = width, height = height, pointsize = pointsize, family = "Courier")
par(mfrow = c(1, 1))
##Obs
i <- as.image.SpatialGridDataFrame(surf.pos1)
plot(t(coords_ord[, S]), typ="n", cex=0.5, xlim=xlim, axes=FALSE, ylab="y", xlab="x") 
     #main = "true")
axis(2, las=1)
axis(1)
image.plot(i, add=TRUE, col=rev(col.pal(length(surf.brks)-1)), zlim=zlim)
dev.off()

png(paste("./pics/map-w2_incp-posm.png", sep = ""), 
    width = width, height = height, pointsize = pointsize, family = "Courier")
par(mfrow = c(1, 1))
##Obs
i <- as.image.SpatialGridDataFrame(surf.pos2)
plot(t(coords_ord[, S]), typ="n", cex=0.5, xlim=xlim, axes=FALSE, ylab="y", xlab="x") 
     #main = "true")
axis(2, las=1)
axis(1)
image.plot(i, add=TRUE, col=rev(col.pal(length(surf.brks)-1)), zlim=zlim)
dev.off()

"""
Out[55]:
RObject{IntSxp}
null device 
          1 

![w1] ![w1_pos] ![w2] ![w2_pos]